linalg_basic Module


Uses


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), public, parameter :: LA_HERMITIAN_TRANSPOSE = 2

Defines a Hermitian transpose operation for a complex-valued matrix.

integer(kind=int32), public, parameter :: LA_NO_OPERATION = 0

Defines no operation should be performed on the matrix.

integer(kind=int32), public, parameter :: LA_TRANSPOSE = 1

Defines a transpose operation.


Interfaces

public interface band_diag_mtx_mult

An interface to the banded diagonal matrix multiplication routines.

  • private subroutine band_diag_mtx_mult_dbl(left, m, kl, ku, alpha, a, b, err)

    Performs the matrix operation or where is a banded matrix and is a diagonal matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: left

    A logical flag indicating whether to perform the operation (TRUE) or (FALSE).

    integer(kind=int32), intent(in) :: m

    The number of rows in the banded matrix .

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    real(kind=real64), intent(inout), dimension(:,:) :: a

    The banded matrix to multiply.

    real(kind=real64), intent(in), dimension(:) :: b

    The diagonal matrix to multiply by.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine band_diag_mtx_mult_cmplx(left, m, kl, ku, alpha, a, b, err)

    Performs the matrix operation or where is a banded matrix and is a diagonal matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: left

    A logical flag indicating whether to perform the operation (TRUE) or (FALSE).

    integer(kind=int32), intent(in) :: m

    The number of rows in the banded matrix .

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    complex(kind=real64), intent(inout), dimension(:,:) :: a

    The banded matrix to multiply.

    complex(kind=real64), intent(in), dimension(:) :: b

    The diagonal matrix to multiply by.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface band_mtx_mult

An interface to the banded matrix multiplication routines.

  • private subroutine band_mtx_vec_mult_dbl(trans, kl, ku, alpha, a, x, beta, y, err)

    Performs the matrix operation or where is a banded matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: trans

    A logical flag indicating whether to perform the operation (FALSE) or (TRUE).

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix .

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix .

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    real(kind=real64), intent(in), dimension(:,:) :: a

    The banded matrix to multiply by.

    real(kind=real64), intent(in), dimension(:) :: x

    The vector to multiply by.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply by.

    real(kind=real64), intent(inout), dimension(:) :: y

    On input, the vector to multiply. On output, the result of the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine band_mtx_vec_mult_cmplx(trans, kl, ku, alpha, a, x, beta, y, err)

    Performs the matrix operation where is a banded matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: trans

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix .

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix .

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    complex(kind=real64), intent(in), dimension(:,:) :: a

    The banded matrix to multiply by.

    complex(kind=real64), intent(in), dimension(:) :: x

    The vector to multiply by.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply by.

    complex(kind=real64), intent(inout), dimension(:) :: y

    On input, the vector to multiply. On output, the result of the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface band_mtx_to_full_mtx

An interface to the banded matrix to full matrix conversion routines.

  • private subroutine band_to_full_mtx_dbl(kl, ku, b, f, err)

    Converts a banded matrix to a full matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix.

    real(kind=real64), intent(in), dimension(:,:) :: b

    The banded matrix to convert.

    real(kind=real64), intent(out), dimension(:,:) :: f

    The full matrix to store the result in.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine band_to_full_mtx_cmplx(kl, ku, b, f, err)

    Converts a banded matrix to a full matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals in the banded matrix.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals in the banded matrix.

    complex(kind=real64), intent(in), dimension(:,:) :: b

    The banded matrix to convert.

    complex(kind=real64), intent(out), dimension(:,:) :: f

    The full matrix to store the result in.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface banded_to_dense

An interface to the banded to dense matrix conversion routines.

  • private subroutine banded_to_dense_dbl(m, kl, ku, a, x, err)

    Converts a banded matrix to a dense matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: m

    The M-by-N dense matrix.

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals. Must be at least 0.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals. Must be at least 0.

    real(kind=real64), intent(in), dimension(:,:) :: a

    The (KL+KU+1)-by-N banded matrix.

    real(kind=real64), intent(out), dimension(:,:) :: x

    The M-by-N dense matrix.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine banded_to_dense_cmplx(m, kl, ku, a, x, err)

    Converts a banded matrix to a dense matrix.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: m

    The M-by-N dense matrix.

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals. Must be at least 0.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals. Must be at least 0.

    complex(kind=real64), intent(in), dimension(:,:) :: a

    The (KL+KU+1)-by-N banded matrix.

    complex(kind=real64), intent(out), dimension(:,:) :: x

    The M-by-N dense matrix.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface dense_to_banded

An interface to the dense to banded matrix conversion routines.

  • private subroutine dense_to_banded_dbl(a, kl, ku, x, err)

    Converts a banded matrix stored in dense format to a compressed form.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in), dimension(:,:) :: a

    The matrix to convert.

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals. Must be at least 0.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals. Must be at least 0.

    real(kind=real64), intent(out), dimension(:,:) :: x

    The (KL+KU+1)-by-N banded matrix.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine dense_to_banded_cmplx(a, kl, ku, x, err)

    Converts a banded matrix stored in dense format to a compressed form.

    The banded matrix is stored in a compressed form supplied column by column. The following code segment transfers between a full matrix to the bonded matrix storage scheme. \code{fortran} do j = 1, n k = ku + 1 - j do i = max(1, j - ku), min(n, j + kl) a(k + i, j) = matrix(i, j) end do end do \endcode

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(in), dimension(:,:) :: a

    The matrix to convert.

    integer(kind=int32), intent(in) :: kl

    The number of subdiagonals. Must be at least 0.

    integer(kind=int32), intent(in) :: ku

    The number of superdiagonals. Must be at least 0.

    complex(kind=real64), intent(out), dimension(:,:) :: x

    The (KL+KU+1)-by-N banded matrix.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface det

An interface to the determinant routines.

  • private function det_dbl(a, iwork, err) result(x)

    Computes the determinant of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout), dimension(:,:) :: a

    On input, the matrix on which to operate. On output, the LU factored matrix in the form [L\U] where L is unit lower triangular and U is upper triangular. The unit diagonal elements of L are not stored.

    integer(kind=int32), intent(out), optional, target, dimension(:) :: iwork

    An MIN(M, N)-element array used to track row-pivot operations. The array stored pivot information such that row I is interchanged with row IPVT(I). If not supplied, this array is allocated within.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

    Return Value real(kind=real64)

    The determinant of the matrix.

  • private function det_cmplx(a, iwork, err) result(x)

    Computes the determinant of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(inout), dimension(:,:) :: a

    On input, the matrix on which to operate. On output, the LU factored matrix in the form [L\U] where L is unit lower triangular and U is upper triangular. The unit diagonal elements of L are not stored.

    integer(kind=int32), intent(out), optional, target, dimension(:) :: iwork

    An MIN(M, N)-element array used to track row-pivot operations. The array stored pivot information such that row I is interchanged with row IPVT(I). If not supplied, this array is allocated within.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

    Return Value complex(kind=real64)

    The determinant of the matrix.

public interface diag_mtx_mult

An interface to the diagonal matrix multiplication routines.

  • private subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    logical, intent(in) :: trans

    A logical flag indicating if the matrix should be transposed.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    real(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    real(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    real(kind=real64), intent(inout), dimension(:,:) :: b

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    logical, intent(in) :: trans

    A logical flag indicating if the matrix should be transposed.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    real(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    complex(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    integer(kind=int32), intent(in) :: opb

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    complex(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    complex(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    integer(kind=int32), intent(in) :: opb

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    complex(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    complex(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    complex(kind=real64), intent(inout), dimension(:,:) :: b

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    integer(kind=int32), intent(in) :: opb

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    complex(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    complex(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)

    Performs the matrix operation or where is a diagonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside

    A logical flag indicating if the diagonal matrix is on the left.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:) :: a

    The diagonal matrix in the operation.

    complex(kind=real64), intent(inout), dimension(:,:) :: b

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine diag_mtx_sparse_mult(lside, alpha, a, b, err)

    Performs the matrix operation or where is a diagonal matrix and is a sparse matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: lside
    real(kind=real64), intent(in) :: alpha
    real(kind=real64), intent(in), dimension(:) :: a
    class(csr_matrix), intent(inout) :: b
    class(errors), intent(inout), optional, target :: err

public interface extract_diagonal

An interface to the diagonal extraction routines.

  • private subroutine extract_diagonal_dbl(a, diag, err)

    Extracts the diagonal of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in), dimension(:,:) :: a

    The M-by-N matrix.

    real(kind=real64), intent(out), dimension(:) :: diag

    The MIN(M, N) element array for the diagonal elements.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine extract_diagonal_cmplx(a, diag, err)

    Extracts the diagonal of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(in), dimension(:,:) :: a

    The M-by-N matrix.

    complex(kind=real64), intent(out), dimension(:) :: diag

    The MIN(M, N) element array for the diagonal elements.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine extract_diagonal_csr(a, diag, err)

    Extracts the diagonal of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    class(csr_matrix), intent(in) :: a

    The M-by-N matrix.

    real(kind=real64), intent(out), dimension(:) :: diag

    The MIN(M, N) element array for the diagonal elements.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface mtx_mult

An interface to the matrix multiplication routines.

  • private subroutine mtx_mult_mtx(transa, transb, alpha, a, b, beta, c, err)

    Performs the matrix operation .

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: transa

    A logical flag indicating if the matrix should be transposed.

    logical, intent(in) :: transb

    A logical flag indicating if the matrix should be transposed.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:,:) :: a

    The matrix in the operation.

    real(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    real(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)

    Performs the matrix-vector operation .

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: trans

    A logical flag indicating if the matrix should be transposed.

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    real(kind=real64), intent(in), dimension(:,:) :: a

    The matrix in the operation.

    real(kind=real64), intent(in), dimension(:) :: b

    The vector in the operation.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply the vector .

    real(kind=real64), intent(inout), dimension(:) :: c

    The vector in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)

    Performs the matrix operation .

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: opa

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    integer(kind=int32), intent(in) :: opb

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:,:) :: a

    The matrix in the operation.

    complex(kind=real64), intent(in), dimension(:,:) :: b

    The matrix in the operation.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply the matrix .

    complex(kind=real64), intent(inout), dimension(:,:) :: c

    The matrix in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)

    Performs the matrix-vector operation .

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: opa

    An integer flag indicating the operation to perform on matrix . Possible options are:

    • LA_NO_OPERATION: No operation is performed on matrix.

    • LA_TRANSPOSE: The transpose of matrix is used.

    • LA_HERMITIAN_TRANSPOSE: The Hermitian transpose of matrix is used.

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the product of and .

    complex(kind=real64), intent(in), dimension(:,:) :: a

    The matrix in the operation.

    complex(kind=real64), intent(in), dimension(:) :: b

    The vector in the operation.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply the vector .

    complex(kind=real64), intent(inout), dimension(:) :: c

    The vector in the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface mtx_rank

An interface to the matrix rank routines.

  • private function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)

    Computes the rank of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout), dimension(:,:) :: a

    The matrix.

    real(kind=real64), intent(in), optional :: tol

    An optional input, that if supplied, overrides the default tolerance on singular values such that singular values less than this tolerance are treated as zero. The default tolerance is: MAX(M, N) * EPS * MAX(S). If the supplied value is less than the smallest value that causes an overflow if inverted, the tolerance reverts back to its default value, and the operation continues; however, a warning message is issued.

    real(kind=real64), intent(out), optional, target, dimension(:) :: work

    An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork. If not provided, the memory required is allocated within.

    integer(kind=int32), intent(out), optional :: olwork

    An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

    Return Value integer(kind=int32)

    The rank of the matrix.

  • private function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)

    Computes the rank of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(inout), dimension(:,:) :: a

    The matrix.

    real(kind=real64), intent(in), optional :: tol

    An optional input, that if supplied, overrides the default tolerance on singular values such that singular values less than this tolerance are treated as zero. The default tolerance is: MAX(M, N) * EPS * MAX(S). If the supplied value is less than the smallest value that causes an overflow if inverted, the tolerance reverts back to its default value, and the operation continues; however, a warning message is issued.

    complex(kind=real64), intent(out), optional, target, dimension(:) :: work

    An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork. If not provided, the memory required is allocated within.

    integer(kind=int32), intent(out), optional :: olwork

    An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

    real(kind=real64), intent(out), optional, target, dimension(:) :: rwork

    An optional input, that if provided, prevents any local memory allocation for real-valued workspace arrays. If not provided, the memory required is allocated within. If provided, the length of the array must be at least 6 * MIN(M, N).

    class(errors), intent(inout), optional, target :: err

    The rank of the matrix.

    Return Value integer(kind=int32)

    The rank of the matrix.

public interface rank1_update

An interface to the rank-1 update routines.

  • private subroutine rank1_update_dbl(alpha, x, y, a, err)

    Performs a rank-1 update of a matrix of the form .

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in) :: alpha

    The scalar to multiply the outer product of and .

    real(kind=real64), intent(in), dimension(:) :: x

    The vector in the outer product.

    real(kind=real64), intent(in), dimension(:) :: y

    The vector in the outer product.

    real(kind=real64), intent(inout), dimension(:,:) :: a

    The matrix to update.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine rank1_update_cmplx(alpha, x, y, a, err)

    Performs a rank-1 update of a matrix of the form .

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply the outer product of and .

    complex(kind=real64), intent(in), dimension(:) :: x

    The vector in the outer product.

    complex(kind=real64), intent(in), dimension(:) :: y

    The vector in the outer product.

    complex(kind=real64), intent(inout), dimension(:,:) :: a

    The matrix to update.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface recip_mult_array

An interface to the reciprocal multiplication routines.

  • private subroutine recip_mult_array_dbl(a, x)

    Computes the product of a scalar and a vector, where the scalar is the reciprocal of the scalar A.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in) :: a

    The scalar A, which is the reciprocal of the scalar to multiply by.

    real(kind=real64), intent(inout), dimension(:) :: x

    On input, the vector to multiply. On output, the product of the vector and the scalar reciprocal.

public interface swap

An interface to the swap routines.

  • private subroutine swap_dbl(x, y, err)

    Swaps the contents of two arrays.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout), dimension(:) :: x

    On input, the first array to swap. On output, the contents of the first array are copied to the second array.

    real(kind=real64), intent(inout), dimension(:) :: y

    On input, the second array to swap. On output, the contents of the second array are copied to the first array.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine swap_cmplx(x, y, err)

    Swaps the contents of two arrays.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(inout), dimension(:) :: x

    On input, the first array to swap. On output, the contents of the first array are copied to the second array.

    complex(kind=real64), intent(inout), dimension(:) :: y

    On input, the second array to swap. On output, the contents of the second array are copied to the first array.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

public interface trace

An interface to the trace routines.

  • private pure function trace_dbl(x) result(y)

    Computes the trace of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in), dimension(:,:) :: x

    The matrix.

    Return Value real(kind=real64)

    The trace of the matrix.

  • private pure function trace_cmplx(x) result(y)

    Computes the trace of a matrix.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(in), dimension(:,:) :: x

    The matrix.

    Return Value complex(kind=real64)

    The trace of the matrix.

public interface tri_mtx_mult

An interface to the triangular matrix multiplication routines.

  • private subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)

    Performs the matrix operation or where is a triangular matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: upper

    A logical flag indicating whether the matrix A is upper triangular (TRUE) or lower triangular (FALSE).

    real(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    real(kind=real64), intent(in), dimension(:,:) :: a

    The triangular matrix to multiply by.

    real(kind=real64), intent(in) :: beta

    The scalar to multiply by.

    real(kind=real64), intent(inout), dimension(:,:) :: b

    On input, the matrix to multiply. On output, the result of the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.

  • private subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)

    Performs the matrix operation or where is a triangular matrix.

    Arguments

    Type IntentOptional Attributes Name
    logical, intent(in) :: upper

    A logical flag indicating whether the matrix A is upper triangular (TRUE) or lower triangular (FALSE).

    complex(kind=real64), intent(in) :: alpha

    The scalar to multiply by.

    complex(kind=real64), intent(in), dimension(:,:) :: a

    The triangular matrix to multiply by.

    complex(kind=real64), intent(in) :: beta

    The scalar to multiply by.

    complex(kind=real64), intent(inout), dimension(:,:) :: b

    On input, the matrix to multiply. On output, the result of the operation.

    class(errors), intent(inout), optional, target :: err

    An error object to report any errors that occur.